Source of the materials: Biopython Tutorial and Cookbook (adapted)

Sequence Input/Output

In this noteboo we’ll discuss in more detail the Bio.SeqIO module, which was briefly introduced before. This aims to provide a simple interface for working with assorted sequence file formats in a uniform way. See also the Bio.SeqIO wiki page (http://biopython.org/wiki/SeqIO), and the built in documentation (also http://biopython.org/DIST/docs/api/Bio.SeqIO-module.html :

In [1]:
import gzip
from io import StringIO

from Bio import Entrez
from Bio import ExPASy
from Bio import SeqIO
from Bio.Alphabet import generic_protein
from Bio.Seq import Seq
from Bio.SeqRecord import SeqRecord
from Bio.SeqUtils.CheckSum import seguid

The catch is that you have to work with SeqRecord objects plus annotation like an identifier and description.

Parsing or Reading Sequences

The workhorse function Bio.SeqIO.parse() is used to read in sequence data as SeqRecord objects. This function expects two arguments:

  • The first argument is a handle to read the data from, or a filename. A handle is typically a file opened for reading, but could be the output from a command line program, or data downloaded from the internet.
  • The second argument is a lower case string specifying sequence format – we don’t try and guess the file format for you! See http://biopython.org/wiki/SeqIO for a full listing of supported formats.

There is an optional argument alphabet to specify the alphabet to be used. This is useful for file formats like FASTA where otherwise Bio.SeqIO will default to a generic alphabet.

The Bio.SeqIO.parse() function returns an iterator which gives SeqRecord objects. Iterators are typically used in a for loop as shown below.

Sometimes you’ll find yourself dealing with files which contain only a single record. For this situation use the function Bio.SeqIO.read() which takes the same arguments. Provided there is one and only one record in the file, this is returned as a SeqRecord object. Otherwise an exception is raised.

Reading Sequence Files

In general Bio.SeqIO.parse() is used to read in sequence files as SeqRecord objects, and is typically used with a for loop like this:

In [2]:
# we show the first 3 only
for i, seq_record in enumerate(SeqIO.parse("data/ls_orchid.fasta", "fasta")):
    print(seq_record.id)
    print(repr(seq_record.seq))
    print(len(seq_record))
    if i == 2:
        break
gi|2765658|emb|Z78533.1|CIZ78533
Seq('CGTAACAAGGTTTCCGTAGGTGAACCTGCGGAAGGATCATTGATGAGACCGTGG...CGC', SingleLetterAlphabet())
740
gi|2765657|emb|Z78532.1|CCZ78532
Seq('CGTAACAAGGTTTCCGTAGGTGAACCTGCGGAAGGATCATTGTTGAGACAACAG...GGC', SingleLetterAlphabet())
753
gi|2765656|emb|Z78531.1|CFZ78531
Seq('CGTAACAAGGTTTCCGTAGGTGAACCTGCGGAAGGATCATTGTTGAGACAGCAG...TAA', SingleLetterAlphabet())
748

The above example is repeated from a previous notebook, and will load the orchid DNA sequences in the FASTA format file. If instead you wanted to load a GenBank format file then all you need to do is change the filename and the format string:

In [3]:
#we show the frist 3
for i, seq_record in enumerate(SeqIO.parse("data/ls_orchid.gbk", "genbank")):
    print(seq_record.id)
    print(seq_record.seq)
    print(len(seq_record))
    if i == 2:
        break
Z78533.1
CGTAACAAGGTTTCCGTAGGTGAACCTGCGGAAGGATCATTGATGAGACCGTGGAATAAACGATCGAGTGAATCCGGAGGACCGGTGTACTCAGCTCACCGGGGGCATTGCTCCCGTGGTGACCCTGATTTGTTGTTGGGCCGCCTCGGGAGCGTCCATGGCGGGTTTGAACCTCTAGCCCGGCGCAGTTTGGGCGCCAAGCCATATGAAAGCATCACCGGCGAATGGCATTGTCTTCCCCAAAACCCGGAGCGGCGGCGTGCTGTCGCGTGCCCAATGAATTTTGATGACTCTCGCAAACGGGAATCTTGGCTCTTTGCATCGGATGGAAGGACGCAGCGAAATGCGATAAGTGGTGTGAATTGCAAGATCCCGTGAACCATCGAGTCTTTTGAACGCAAGTTGCGCCCGAGGCCATCAGGCTAAGGGCACGCCTGCTTGGGCGTCGCGCTTCGTCTCTCTCCTGCCAATGCTTGCCCGGCATACAGCCAGGCCGGCGTGGTGCGGATGTGAAAGATTGGCCCCTTGTGCCTAGGTGCGGCGGGTCCAAGAGCTGGTGTTTTGATGGCCCGGAACCCGGCAAGAGGTGGACGGATGCTGGCAGCAGCTGCCGTGCGAATCCCCCATGTTGTCGTGCTTGTCGGACAGGCAGGAGAACCCTTCCGAACCCCAATGGAGGGCGGTTGACCGCCATTCGGATGTGACCCCAGGTCAGGCGGGGGCACCCGCTGAGTTTACGC
740
Z78532.1
CGTAACAAGGTTTCCGTAGGTGAACCTGCGGAAGGATCATTGTTGAGACAACAGAATATATGATCGAGTGAATCTGGAGGACCTGTGGTAACTCAGCTCGTCGTGGCACTGCTTTTGTCGTGACCCTGCTTTGTTGTTGGGCCTCCTCAAGAGCTTTCATGGCAGGTTTGAACTTTAGTACGGTGCAGTTTGCGCCAAGTCATATAAAGCATCACTGATGAATGACATTATTGTCAGAAAAAATCAGAGGGGCAGTATGCTACTGAGCATGCCAGTGAATTTTTATGACTCTCGCAACGGATATCTTGGCTCTAACATCGATGAAGAACGCAGCTAAATGCGATAAGTGGTGTGAATTGCAGAATCCCGTGAACCATCGAGTCTTTGAACGCAAGTTGCGCTCGAGGCCATCAGGCTAAGGGCACGCCTGCCTGGGCGTCGTGTGTTGCGTCTCTCCTACCAATGCTTGCTTGGCATATCGCTAAGCTGGCATTATACGGATGTGAATGATTGGCCCCTTGTGCCTAGGTGCGGTGGGTCTAAGGATTGTTGCTTTGATGGGTAGGAATGTGGCACGAGGTGGAGAATGCTAACAGTCATAAGGCTGCTATTTGAATCCCCCATGTTGTTGTATTTTTTCGAACCTACACAAGAACCTAATTGAACCCCAATGGAGCTAAAATAACCATTGGGCAGTTGATTTCCATTCAGATGCGACCCCAGGTCAGGCGGGGCCACCCGCTGAGTTGAGGC
753
Z78531.1
CGTAACAAGGTTTCCGTAGGTGAACCTGCGGAAGGATCATTGTTGAGACAGCAGAACATACGATCGAGTGAATCCGGAGGACCCGTGGTTACACGGCTCACCGTGGCTTTGCTCTCGTGGTGAACCCGGTTTGCGACCGGGCCGCCTCGGGAACTTTCATGGCGGGTTTGAACGTCTAGCGCGGCGCAGTTTGCGCCAAGTCATATGGAGCGTCACCGATGGATGGCATTTTTGTCAAGAAAAACTCGGAGGGGCGGCGTCTGTTGCGCGTGCCAATGAATTTATGACGACTCTCGGCAACGGGATATCTGGCTCTTGCATCGATGAAGAACGCAGCGAAATGCGATAAGTGGTGTGAATTGCAGAATCCCGCGAACCATCGAGTCTTTGAACGCAAGTTGCGCCCGAGGCCATCAGGCTAAGGGCACGCCTGCCTGGGCGTCGTGTGCTGCGTCTCTCCTGATAATGCTTGATTGGCATGCGGCTAGTCTGTCATTGTGAGGACGTGAAAGATTGGCCCCTTGCGCCTAGGTGCGGCGGGTCTAAGCATCGGTGTTCTGATGGCCCGGAACTTGGCAGTAGGTGGAGGATGCTGGCAGCCGCAAGGCTGCCGTTCGAATCCCCCGTGTTGTCGTACTCGTCAGGCCTACAGAAGAACCTGTTTGAACCCCCAGTGGACGCAAAACCGCCCTCGGGCGGTGATTTCCATTCAGATGCGACCCCAGTCAGGCGGGCCACCCGTGAGTAA
748

Similarly, if you wanted to read in a file in another file format, then assuming Bio.SeqIO.parse() supports it you would just need to change the format string as appropriate, for example ‘swiss’ for SwissProt files or ‘embl’ for EMBL text files. There is a full listing on the wiki page (http://biopython.org/wiki/SeqIO) and in the built in documentation (also :raw-latex:`\http`://biopython.org/DIST/docs/api/Bio.SeqIO-module.html).

Another very common way to use a Python iterator is within a list comprehension (or a generator expression). For example, if all you wanted to extract from the file was a list of the record identifiers we can easily do this with the following list comprehension:

In [4]:
[seq_record.id for seq_record in SeqIO.parse("data/ls_orchid.gbk", "genbank")][:10]  # ten only
Out[4]:
['Z78533.1',
 'Z78532.1',
 'Z78531.1',
 'Z78530.1',
 'Z78529.1',
 'Z78527.1',
 'Z78526.1',
 'Z78525.1',
 'Z78524.1',
 'Z78523.1']

Iterating over the records in a sequence file

In the above examples, we have usually used a for loop to iterate over all the records one by one. You can use the for loop with all sorts of Python objects (including lists, tuples and strings) which support the iteration interface.

The object returned by Bio.SeqIO is actually an iterator which returns SeqRecord objects. You get to see each record in turn, but once and only once. The plus point is that an iterator can save you memory when dealing with large files.

Instead of using a for loop, can also use the next() function on an iterator to step through the entries, like this:

In [5]:
record_iterator = SeqIO.parse("data/ls_orchid.fasta", "fasta")

first_record = next(record_iterator)
print(first_record.id)
print(first_record.description)

second_record = next(record_iterator)
print(second_record.id)
print(second_record.description)
gi|2765658|emb|Z78533.1|CIZ78533
gi|2765658|emb|Z78533.1|CIZ78533 C.irapeanum 5.8S rRNA gene and ITS1 and ITS2 DNA
gi|2765657|emb|Z78532.1|CCZ78532
gi|2765657|emb|Z78532.1|CCZ78532 C.californicum 5.8S rRNA gene and ITS1 and ITS2 DNA

Note that if you try to use next() and there are no more results, you’ll get the special StopIteration exception.

One special case to consider is when your sequence files have multiple records, but you only want the first one. In this situation the following code is very concise:

In [6]:
next(SeqIO.parse("data/ls_orchid.gbk", "genbank"))
Out[6]:
SeqRecord(seq=Seq('CGTAACAAGGTTTCCGTAGGTGAACCTGCGGAAGGATCATTGATGAGACCGTGG...CGC', IUPACAmbiguousDNA()), id='Z78533.1', name='Z78533', description='C.irapeanum 5.8S rRNA gene and ITS1 and ITS2 DNA.', dbxrefs=[])

A word of warning here – using the next() function like this will silently ignore any additional records in the file. If your files have one and only one record, like some of the online examples later in this chapter, or a GenBank file for a single chromosome, then use the new Bio.SeqIO.read() function instead. This will check there are no extra unexpected records present.

Getting a list of the records in a sequence file

In the previous section we talked about the fact that Bio.SeqIO.parse() gives you a SeqRecord iterator, and that you get the records one by one. Very often you need to be able to access the records in any order. The Python list data type is perfect for this, and we can turn the record iterator into a list of SeqRecord objects using the built-in Python function list() like so:

In [7]:
records = list(SeqIO.parse("data/ls_orchid.gbk", "genbank"))

print("Found %i records" % len(records))

print("The last record")
last_record = records[-1] #using Python's list tricks
print(last_record.id)
print(repr(last_record.seq))
print(len(last_record))

print("The first record")
first_record = records[0] #remember, Python counts from zero
print(first_record.id)
print(repr(first_record.seq))
print(len(first_record))

Found 94 records
The last record
Z78439.1
Seq('CATTGTTGAGATCACATAATAATTGATCGAGTTAATCTGGAGGATCTGTTTACT...GCC', IUPACAmbiguousDNA())
592
The first record
Z78533.1
Seq('CGTAACAAGGTTTCCGTAGGTGAACCTGCGGAAGGATCATTGATGAGACCGTGG...CGC', IUPACAmbiguousDNA())
740

You can of course still use a for loop with a list of SeqRecord objects. Using a list is much more flexible than an iterator (for example, you can determine the number of records from the length of the list), but does need more memory because it will hold all the records in memory at once.

Extracting data

The SeqRecord object and its annotation structures are described more fully in in another notebook. As an example of how annotations are stored, we’ll look at the output from parsing the first record in the orchid GenBank file.

In [8]:
record_iterator = SeqIO.parse("data/ls_orchid.gbk", "genbank")
first_record = next(record_iterator)
print(first_record)
ID: Z78533.1
Name: Z78533
Description: C.irapeanum 5.8S rRNA gene and ITS1 and ITS2 DNA.
Number of features: 5
/organism=Cypripedium irapeanum
/taxonomy=['Eukaryota', 'Viridiplantae', 'Streptophyta', 'Embryophyta', 'Tracheophyta', 'Spermatophyta', 'Magnoliophyta', 'Liliopsida', 'Asparagales', 'Orchidaceae', 'Cypripedioideae', 'Cypripedium']
/sequence_version=1
/keywords=['5.8S ribosomal RNA', '5.8S rRNA gene', 'internal transcribed spacer', 'ITS1', 'ITS2']
/date=30-NOV-2006
/accessions=['Z78533']
/source=Cypripedium irapeanum
/references=[Reference(title='Phylogenetics of the slipper orchids (Cypripedioideae: Orchidaceae): nuclear rDNA ITS sequences', ...), Reference(title='Direct Submission', ...)]
/gi=2765658
/data_file_division=PLN
Seq('CGTAACAAGGTTTCCGTAGGTGAACCTGCGGAAGGATCATTGATGAGACCGTGG...CGC', IUPACAmbiguousDNA())

This gives a human readable summary of most of the annotation data for the SeqRecord. For this example we’re going to use the .annotations attribute which is just a Python dictionary. The contents of this annotations dictionary were shown when we printed the record above. You can also print them out directly:

print(first_record.annotations)

Like any Python dictionary, you can easily get a list of the keys:

print(first_record.annotations.keys())

or values:

print(first_record.annotations.values())

In general, the annotation values are strings, or lists of strings. One special case is any references in the file get stored as reference objects.

Suppose you wanted to extract a list of the species from the GenBank file. The information we want, Cypripedium irapeanum, is held in the annotations dictionary under ‘source’ and ‘organism’, which we can access like this:

In [9]:
print(first_record.annotations["source"])
print(first_record.annotations["organism"])
Cypripedium irapeanum
Cypripedium irapeanum

In general, ‘organism’ is used for the scientific name (in Latin, e.g. Arabidopsis thaliana), while ‘source’ will often be the common name (e.g. thale cress). In this example, as is often the case, the two fields are identical.

Now let’s go through all the records, building up a list of the species each orchid sequence is from:

In [10]:
all_species = []
for seq_record in SeqIO.parse("data/ls_orchid.gbk", "genbank"):
    all_species.append(seq_record.annotations["organism"])
print(all_species[:10])  # we print only 10
['Cypripedium irapeanum', 'Cypripedium californicum', 'Cypripedium fasciculatum', 'Cypripedium margaritaceum', 'Cypripedium lichiangense', 'Cypripedium yatabeanum', 'Cypripedium guttatum', 'Cypripedium acaule', 'Cypripedium formosanum', 'Cypripedium himalaicum']

Another way of writing this code is to use a list comprehension:

In [11]:
all_species = [seq_record.annotations["organism"] for seq_record in \
               SeqIO.parse("data/ls_orchid.gbk", "genbank")]
print(all_species[:10])
['Cypripedium irapeanum', 'Cypripedium californicum', 'Cypripedium fasciculatum', 'Cypripedium margaritaceum', 'Cypripedium lichiangense', 'Cypripedium yatabeanum', 'Cypripedium guttatum', 'Cypripedium acaule', 'Cypripedium formosanum', 'Cypripedium himalaicum']

Great. That was pretty easy because GenBank files are annotated in a standardised way.

Now, let’s suppose you wanted to extract a list of the species from a FASTA file, rather than the GenBank file. The bad news is you will have to write some code to extract the data you want from the record’s description line - if the information is in the file in the first place! Our example FASTA format file starts like this:

>gi|2765658|emb|Z78533.1|CIZ78533 C.irapeanum 5.8S rRNA gene and ITS1 and ITS2 DNA
CGTAACAAGGTTTCCGTAGGTGAACCTGCGGAAGGATCATTGATGAGACCGTGGAATAAACGATCGAGTG
AATCCGGAGGACCGGTGTACTCAGCTCACCGGGGGCATTGCTCCCGTGGTGACCCTGATTTGTTGTTGGG
...

You can check by hand, but for every record the species name is in the description line as the second word. This means if we break up each record’s .description at the spaces, then the species is there as field number one (field zero is the record identifier). That means we can do this:

In [12]:
all_species = []
for seq_record in SeqIO.parse("data/ls_orchid.fasta", "fasta"):
    all_species.append(seq_record.description.split()[1])
print(all_species[:10])
['C.irapeanum', 'C.californicum', 'C.fasciculatum', 'C.margaritaceum', 'C.lichiangense', 'C.yatabeanum', 'C.guttatum', 'C.acaule', 'C.formosanum', 'C.himalaicum']

The concise alternative using list comprehensions would be:

In [13]:
all_species == [seq_record.description.split()[1] for seq_record in \
                SeqIO.parse("data/ls_orchid.fasta", "fasta")]
print(all_species[:10])
['C.irapeanum', 'C.californicum', 'C.fasciculatum', 'C.margaritaceum', 'C.lichiangense', 'C.yatabeanum', 'C.guttatum', 'C.acaule', 'C.formosanum', 'C.himalaicum']

In general, extracting information from the FASTA description line is not very nice. If you can get your sequences in a well annotated file format like GenBank or EMBL, then this sort of annotation information is much easier to deal with.

Parsing sequences from compressed files

In the previous section, we looked at parsing sequence data from a file. Instead of using a filename, you can give Bio.SeqIO a handle, and in this section we’ll use handles to parse sequence from compressed files.

As you’ll have seen above, we can use Bio.SeqIO.read() or Bio.SeqIO.parse() with a filename - for instance this quick example calculates the total length of the sequences in a multiple record GenBank file using a generator expression:

In [14]:
print(sum(len(r) for r in SeqIO.parse("data/ls_orchid.gbk", "gb")))
67518

Here we use a file handle instead, using the :raw-latex:`\verb|with|` statement to close the handle automatically:

In [15]:
with open("data/ls_orchid.gbk") as handle:
    print(sum(len(r) for r in SeqIO.parse(handle, "gb")))
67518

Or, the old fashioned way where you manually close the handle:

In [16]:
handle = open("data/ls_orchid.gbk")
print(sum(len(r) for r in SeqIO.parse(handle, "gb")))
handle.close()
67518

Now, suppose we have a gzip compressed file instead? These are very commonly used on Linux. We can use Python’s gzip module to open the compressed file for reading - which gives us a handle object:

In [17]:
handle = gzip.open("data/ls_orchid.gbk.gz", "rt")
print(sum(len(r) for r in SeqIO.parse(handle, "gb")))
handle.close()
67518

There is a gzip (GNU Zip) variant called BGZF (Blocked GNU Zip Format), which can be treated like an ordinary gzip file for reading, but has advantages for random access later which we’ll talk about later.

Parsing sequences from the net

In the previous sections, we looked at parsing sequence data from a file (using a filename or handle), and from compressed files (using a handle). Here we’ll use Bio.SeqIO with another type of handle, a network connection, to download and parse sequences from the internet.

Note that just because you can download sequence data and parse it into a SeqRecord object in one go doesn’t mean this is a good idea. In general, you should probably download sequences once and save them to a file for reuse.

Parsing GenBank records from the net

Let’s just connect to the NCBI and get a few Opuntia (prickly-pear) sequences from GenBank using their GI numbers.

First of all, let’s fetch just one record. If you don’t care about the annotations and features downloading a FASTA file is a good choice as these are compact. Now remember, when you expect the handle to contain one and only one record, use the Bio.SeqIO.read() function:

In [18]:
Entrez.email = "A.N.Other@example.com"
handle = Entrez.efetch(db="nucleotide", rettype="fasta", retmode="text", id="6273291")
seq_record = SeqIO.read(handle, "fasta")
handle.close()
print("%s with %i features" % (seq_record.id, len(seq_record.features)))
gi|6273291|gb|AF191665.1|AF191665 with 0 features

The NCBI will also let you ask for the file in other formats, in particular as a GenBank file. Until Easter 2009, the Entrez EFetch API let you use ``genbank’’ as the return type, however the NCBI now insist on using the official return types of ‘gb’ (or ‘gp’ for proteins) as described on http://www.ncbi.nlm.nih.gov/entrez/query/static/efetchseq_help.html.

As a result we support ‘gb’ as an alias for ‘genbank’ in Bio.SeqIO.

In [19]:
Entrez.email = "A.N.Other@example.com"
handle = Entrez.efetch(db="nucleotide", rettype="gb", retmode="text", id="6273291")
seq_record = SeqIO.read(handle, "gb") #using "gb" as an alias for "genbank"
handle.close()
print("%s with %i features" % (seq_record.id, len(seq_record.features)))
AF191665.1 with 3 features

Notice this time we have three features.

Now let’s fetch several records. This time the handle contains multiple records, so we must use the :raw-latex:`\verb|Bio.SeqIO.parse()|` function:

In [20]:
Entrez.email = "A.N.Other@example.com"
handle = Entrez.efetch(db="nucleotide", rettype="gb", retmode="text",
                       id="6273291,6273290,6273289")
for seq_record in SeqIO.parse(handle, "gb"):
    print (seq_record.id, seq_record.description[:50] + "...")
    print ("Sequence length %i," % len(seq_record))
    print ("%i features," % len(seq_record.features))
    print ("from: %s" % seq_record.annotations["source"])
handle.close()
AF191665.1 Opuntia marenae rpl16 gene; chloroplast gene for c...
Sequence length 902,
3 features,
from: chloroplast Grusonia marenae
AF191664.1 Opuntia clavata rpl16 gene; chloroplast gene for c...
Sequence length 899,
3 features,
from: chloroplast Grusonia clavata
AF191663.1 Opuntia bradtiana rpl16 gene; chloroplast gene for...
Sequence length 899,
3 features,
from: chloroplast Grusonia bradtiana

Parsing SwissProt sequences from the net

Now let’s use a handle to download a SwissProt file from ExPASy (SwissProt will be covered in detail in another notebook). As mentioned above, when you expect the handle to contain one and only one record, use the Bio.SeqIO.read() function:

In [21]:
handle = ExPASy.get_sprot_raw("O23729")
seq_record = SeqIO.read(handle, "swiss")
handle.close()
print(seq_record.id)
print(seq_record.name)
print(seq_record.description)
print(repr(seq_record.seq))
print("Length %i" % len(seq_record))
print(seq_record.annotations["keywords"])
O23729
CHS3_BROFI
RecName: Full=Chalcone synthase 3; EC=2.3.1.74; AltName: Full=Naringenin-chalcone synthase 3;
Seq('MAPAMEEIRQAQRAEGPAAVLAIGTSTPPNALYQADYPDYYFRITKSEHLTELK...GAE', ProteinAlphabet())
Length 394
['Acyltransferase', 'Flavonoid biosynthesis', 'Transferase']

Sequence files as dictionaries

We’re now going to introduce three related functions in the :raw-latex:`\verb|Bio.SeqIO|` module which allow dictionary like random access to a multi-sequence file. There is a trade off here between flexibility and memory usage. In summary:

  • Bio.SeqIO.to_dict() is the most flexible but also the most memory demanding option. This is basically a helper function to build a normal Python dictionary with each entry held as a SeqRecord object in memory, allowing you to modify the records.
  • Bio.SeqIO.index() is a useful middle ground, acting like a read only dictionary and parsing sequences into SeqRecord objects on demand.
  • Bio.SeqIO.index_db() also acts like a read only dictionary but stores the identifiers and file offsets in a file on disk (as an SQLite3 database), meaning it has very low memory requirements, but will be a little bit slower.

Sequence files as Dictionaries – In memory

The next thing that we’ll do with our ubiquitous orchid files is to show how to index them and access them like a database using the Python dictionary data type. This is very useful for moderately large files where you only need to access certain elements of the file, and makes for a nice quick ‘n dirty database. For dealing with larger files where memory becomes a problem, see below.

You can use the function Bio.SeqIO.to_dict() to make a SeqRecord dictionary (in memory). By default this will use each record’s identifier (i.e. the .id attribute) as the key. Let’s try this using our GenBank file:

In [22]:
orchid_dict = SeqIO.to_dict(SeqIO.parse("data/ls_orchid.gbk", "genbank"))

There is just one required argument for Bio.SeqIO.to_dict(), a list or generator giving SeqRecord objects. Here we have just used the output from the SeqIO.parse function. As the name suggests, this returns a Python dictionary.

Since this variable orchid_dict is an ordinary Python dictionary, we can look at all of the keys we have available:

In [23]:
print(len(orchid_dict))
print(orchid_dict.keys())
94
dict_keys(['Z78459.1', 'Z78496.1', 'Z78501.1', 'Z78443.1', 'Z78514.1', 'Z78452.1', 'Z78442.1', 'Z78439.1', 'Z78526.1', 'Z78527.1', 'Z78532.1', 'Z78471.1', 'Z78444.1', 'Z78511.1', 'Z78489.1', 'Z78502.1', 'Z78480.1', 'Z78453.1', 'Z78505.1', 'Z78475.1', 'Z78445.1', 'Z78507.1', 'Z78515.1', 'Z78446.1', 'Z78462.1', 'Z78461.1', 'Z78476.1', 'Z78523.1', 'Z78521.1', 'Z78440.1', 'Z78483.1', 'Z78449.1', 'Z78469.1', 'Z78458.1', 'Z78506.1', 'Z78470.1', 'Z78510.1', 'Z78477.1', 'Z78499.1', 'Z78519.1', 'Z78495.1', 'Z78474.1', 'Z78491.1', 'Z78487.1', 'Z78479.1', 'Z78450.1', 'Z78513.1', 'Z78492.1', 'Z78456.1', 'Z78467.1', 'Z78504.1', 'Z78531.1', 'Z78481.1', 'Z78530.1', 'Z78447.1', 'Z78460.1', 'Z78503.1', 'Z78500.1', 'Z78464.1', 'Z78498.1', 'Z78472.1', 'Z78454.1', 'Z78465.1', 'Z78473.1', 'Z78482.1', 'Z78490.1', 'Z78497.1', 'Z78508.1', 'Z78463.1', 'Z78466.1', 'Z78524.1', 'Z78518.1', 'Z78520.1', 'Z78509.1', 'Z78455.1', 'Z78516.1', 'Z78451.1', 'Z78485.1', 'Z78486.1', 'Z78512.1', 'Z78468.1', 'Z78488.1', 'Z78457.1', 'Z78529.1', 'Z78494.1', 'Z78522.1', 'Z78525.1', 'Z78484.1', 'Z78533.1', 'Z78441.1', 'Z78448.1', 'Z78493.1', 'Z78478.1', 'Z78517.1'])

If you really want to, you can even look at all the records at once:

In [24]:
list(orchid_dict.values())[:5]  # Ok not all at once...
Out[24]:
[SeqRecord(seq=Seq('CGTAACAAGGTTTCCGTAGGTGAACCTGCGGAAGGATCATTGTTGAGATCACAT...TTT', IUPACAmbiguousDNA()), id='Z78459.1', name='Z78459', description='P.dayanum 5.8S rRNA gene and ITS1 and ITS2 DNA.', dbxrefs=[]),
 SeqRecord(seq=Seq('CGTAACAAGGTTTCCGTAGGTGAACCTGCGGAAGGATCATTGTTGAGATCGCAT...AGC', IUPACAmbiguousDNA()), id='Z78496.1', name='Z78496', description='P.armeniacum 5.8S rRNA gene and ITS1 and ITS2 DNA.', dbxrefs=[]),
 SeqRecord(seq=Seq('CGTAACAAGGTTTCCGTAGGTGAACCTGCGGAAGGATCATTGTTGAGACCGCAA...AGA', IUPACAmbiguousDNA()), id='Z78501.1', name='Z78501', description='P.caudatum 5.8S rRNA gene and ITS1 and ITS2 DNA.', dbxrefs=[]),
 SeqRecord(seq=Seq('CGTAACAAGGTTTCCGTAGGTGAACCTGCGGAAGGATCATTGTTGAGATCACAT...AGG', IUPACAmbiguousDNA()), id='Z78443.1', name='Z78443', description='P.lawrenceanum 5.8S rRNA gene and ITS1 and ITS2 DNA.', dbxrefs=[]),
 SeqRecord(seq=Seq('CGTAACAAGGTTTCCGTAGGTGGACCTTCGGGAGGATCATTTTTGAAGCCCCCA...CTA', IUPACAmbiguousDNA()), id='Z78514.1', name='Z78514', description='P.schlimii 5.8S rRNA gene and ITS1 and ITS2 DNA.', dbxrefs=[])]

We can access a single :raw-latex:`\verb|SeqRecord|` object via the keys and manipulate the object as normal:

In [25]:
seq_record = orchid_dict["Z78475.1"]
print(seq_record.description)
print(repr(seq_record.seq))
P.supardii 5.8S rRNA gene and ITS1 and ITS2 DNA.
Seq('CGTAACAAGGTTTCCGTAGGTGAACCTGCGGAAGGATCATTGTTGAGATCACAT...GGT', IUPACAmbiguousDNA())

So, it is very easy to create an in memory ‘database’ of our GenBank records. Next we’ll try this for the FASTA file instead.

Note that those of you with prior Python experience should all be able to construct a dictionary like this ‘by hand’. However, typical dictionary construction methods will not deal with the case of repeated keys very nicely. Using the Bio.SeqIO.to_dict() will explicitly check for duplicate keys, and raise an exception if any are found.

Specifying the dictionary keys

Using the same code as above, but for the FASTA file instead:

In [26]:
orchid_dict = SeqIO.to_dict(SeqIO.parse("data/ls_orchid.fasta", "fasta"))
print(list(orchid_dict.keys())[:5])
['gi|2765611|emb|Z78486.1|PBZ78486', 'gi|2765630|emb|Z78505.1|PSZ78505', 'gi|2765584|emb|Z78459.1|PDZ78459', 'gi|2765632|emb|Z78507.1|PLZ78507', 'gi|2765644|emb|Z78519.1|CPZ78519']

You should recognise these strings from when we parsed the FASTA file earlier. Suppose you would rather have something else as the keys - like the accession numbers. This brings us nicely to SeqIO.to_dict()’s optional argument key_function, which lets you define what to use as the dictionary key for your records.

First you must write your own function to return the key you want (as a string) when given a SeqRecord object. In general, the details of function will depend on the sort of input records you are dealing with. But for our orchids, we can just split up the record’s identifier using the ‘pipe’ character (the vertical line) and return the fourth entry (field three):

In [27]:
def get_accession(record):
    """"Given a SeqRecord, return the accession number as a string.

    e.g. "gi|2765613|emb|Z78488.1|PTZ78488" -> "Z78488.1"
    """
    parts = record.id.split("|")
    assert len(parts) == 5 and parts[0] == "gi" and parts[2] == "emb"
    return parts[3]

Then we can give this function to the SeqIO.to_dict() function to use in building the dictionary:

In [28]:
orchid_dict = SeqIO.to_dict(SeqIO.parse("data/ls_orchid.fasta", "fasta"), key_function=get_accession)
print(orchid_dict.keys())
dict_keys(['Z78459.1', 'Z78496.1', 'Z78501.1', 'Z78443.1', 'Z78514.1', 'Z78452.1', 'Z78442.1', 'Z78439.1', 'Z78526.1', 'Z78527.1', 'Z78532.1', 'Z78471.1', 'Z78444.1', 'Z78511.1', 'Z78489.1', 'Z78502.1', 'Z78480.1', 'Z78453.1', 'Z78505.1', 'Z78475.1', 'Z78445.1', 'Z78507.1', 'Z78515.1', 'Z78446.1', 'Z78462.1', 'Z78461.1', 'Z78476.1', 'Z78523.1', 'Z78521.1', 'Z78440.1', 'Z78483.1', 'Z78449.1', 'Z78469.1', 'Z78458.1', 'Z78506.1', 'Z78470.1', 'Z78510.1', 'Z78477.1', 'Z78499.1', 'Z78519.1', 'Z78495.1', 'Z78474.1', 'Z78491.1', 'Z78487.1', 'Z78479.1', 'Z78450.1', 'Z78513.1', 'Z78492.1', 'Z78456.1', 'Z78467.1', 'Z78504.1', 'Z78531.1', 'Z78481.1', 'Z78530.1', 'Z78447.1', 'Z78460.1', 'Z78503.1', 'Z78500.1', 'Z78464.1', 'Z78498.1', 'Z78472.1', 'Z78454.1', 'Z78465.1', 'Z78473.1', 'Z78482.1', 'Z78490.1', 'Z78497.1', 'Z78508.1', 'Z78463.1', 'Z78466.1', 'Z78524.1', 'Z78518.1', 'Z78520.1', 'Z78509.1', 'Z78455.1', 'Z78516.1', 'Z78451.1', 'Z78485.1', 'Z78486.1', 'Z78512.1', 'Z78468.1', 'Z78488.1', 'Z78457.1', 'Z78529.1', 'Z78494.1', 'Z78522.1', 'Z78525.1', 'Z78484.1', 'Z78533.1', 'Z78441.1', 'Z78448.1', 'Z78493.1', 'Z78478.1', 'Z78517.1'])

Finally, as desired, the new dictionary keys:

Indexing a dictionary using the SEGUID checksum

To give another example of working with dictionaries of SeqRecord objects, we’ll use the SEGUID checksum function. This is a relatively recent checksum, and collisions should be very rare (i.e. two different sequences with the same checksum), an improvement on the CRC64 checksum.

Once again, working with the orchids GenBank file:

In [29]:
for i, record in enumerate(SeqIO.parse("data/ls_orchid.gbk", "genbank")):
    print(record.id, seguid(record.seq))
    if i == 4:  # OK, 5 is enough!
        break
Z78533.1 JUEoWn6DPhgZ9nAyowsgtoD9TTo
Z78532.1 MN/s0q9zDoCVEEc+k/IFwCNF2pY
Z78531.1 xN45pACrTnmBH8a8Y9cWSgoLrwE
Z78530.1 yMhI5UUQfFOPcoJXb9B19XUyYlY
Z78529.1 s1Pnjq9zoSHoI/CG9jQr4GyeMZY

Now, recall the Bio.SeqIO.to_dict() function’s key_function argument expects a function which turns a SeqRecord into a string. We can’t use the seguid() function directly because it expects to be given a Seq object (or a string). However, we can use Python’s lambda feature to create a ‘one off’ function to give to Bio.SeqIO.to_dict() instead:

In [30]:
seguid_dict = SeqIO.to_dict(SeqIO.parse("data/ls_orchid.gbk", "genbank"),
                            lambda rec : seguid(rec.seq))
record = seguid_dict["MN/s0q9zDoCVEEc+k/IFwCNF2pY"]
print(record.id)
print(record.description)
Z78532.1
C.californicum 5.8S rRNA gene and ITS1 and ITS2 DNA.

That should have retrieved the record Z78532.1, the second entry in the file.

Sequence files as Dictionaries – Indexed files

As the previous couple of examples tried to illustrate, using Bio.SeqIO.to_dict() is very flexible. However, because it holds everything in memory, the size of file you can work with is limited by your computer’s RAM. In general, this will only work on small to medium files.

For larger files you should consider Bio.SeqIO.index(), which works a little differently. Although it still returns a dictionary like object, this does not keep everything in memory. Instead, it just records where each record is within the file – when you ask for a particular record, it then parses it on demand.

As an example, let’s use the same GenBank file as before:

In [31]:
orchid_dict = SeqIO.index("data/ls_orchid.gbk", "genbank")
print(len(orchid_dict))
print(orchid_dict.keys())
seq_record = orchid_dict["Z78475.1"]
print(seq_record.description)
print(seq_record.seq)
94
<dict_keyiterator object at 0x7f92241b0db8>
P.supardii 5.8S rRNA gene and ITS1 and ITS2 DNA.
CGTAACAAGGTTTCCGTAGGTGAACCTGCGGAAGGATCATTGTTGAGATCACATAATAATTGATCGAGTTAATCTGGAGGATCAGTTTACTTTGGTCACCCATGGGCATCTGCTCTTGCAGTGACCTGGATTTGCCATCGAGCCTCCTTGGGAGCTTTCTTGCTGGCGATCTAAACCCGTCCCGGCGCAGTTTTGCGCCAAGTCATATGACACATAATTGGAAGGGGGTGGCATGCTGCCTTGACCCTCCCCAAATTATTTTTTTGACAACTCTCAGCAACGGATATCTCGGCTCTTNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNNATCAGGCCAAGGGCACGCCTGCCTGGGCATTGCGAGTCATATCTCTCCCTTAATGAGGCTGTCCATACATACTGTTCAGCCAATGCGGATGTGAGTTTGGCCCCTTGTTCTTTGGTACGGGGGGTCTAAGAGCTGCATGGGCTTTTGATGGTCCAAAATACGGCAAGAGGTGGACGAACTATGCTACAACAAAATTGTTGTGCGAATGCCCCGGGTTGTCGTATTAGATGGGCCAGCATAATCTAAAGACCCTTTTGAACCCCATTGGAGGCCCATCAACCCATGATCAGTTGACGGCCATTTGGTTGCGACCCAGGTCAGGT

Note that Bio.SeqIO.index() won’t take a handle, but only a filename. There are good reasons for this, but it is a little technical. The second argument is the file format (a lower case string as used in the other Bio.SeqIO functions). You can use many other simple file formats, including FASTA and FASTQ files. However, alignment formats like PHYLIP or Clustal are not supported. Finally as an optional argument you can supply an alphabet, or a key function.

Here is the same example using the FASTA file - all we change is the filename and the format name:

In [32]:
orchid_dict = SeqIO.index("data/ls_orchid.fasta", "fasta")
print(len(orchid_dict))
print(list(orchid_dict.keys()))
94
['gi|2765611|emb|Z78486.1|PBZ78486', 'gi|2765630|emb|Z78505.1|PSZ78505', 'gi|2765584|emb|Z78459.1|PDZ78459', 'gi|2765632|emb|Z78507.1|PLZ78507', 'gi|2765644|emb|Z78519.1|CPZ78519', 'gi|2765645|emb|Z78520.1|CSZ78520', 'gi|2765650|emb|Z78525.1|CAZ78525', 'gi|2765571|emb|Z78446.1|PAZ78446', 'gi|2765575|emb|Z78450.1|PPZ78450', 'gi|2765615|emb|Z78490.1|PFZ78490', 'gi|2765577|emb|Z78452.1|PBZ78452', 'gi|2765656|emb|Z78531.1|CFZ78531', 'gi|2765601|emb|Z78476.1|PGZ78476', 'gi|2765588|emb|Z78463.1|PGZ78463', 'gi|2765633|emb|Z78508.1|PLZ78508', 'gi|2765600|emb|Z78475.1|PSZ78475', 'gi|2765634|emb|Z78509.1|PPZ78509', 'gi|2765655|emb|Z78530.1|CMZ78530', 'gi|2765616|emb|Z78491.1|PCZ78491', 'gi|2765631|emb|Z78506.1|PLZ78506', 'gi|2765620|emb|Z78495.1|PEZ78495', 'gi|2765579|emb|Z78454.1|PFZ78454', 'gi|2765566|emb|Z78441.1|PSZ78441', 'gi|2765581|emb|Z78456.1|PTZ78456', 'gi|2765612|emb|Z78487.1|PHZ78487', 'gi|2765626|emb|Z78501.1|PCZ78501', 'gi|2765582|emb|Z78457.1|PCZ78457', 'gi|2765568|emb|Z78443.1|PLZ78443', 'gi|2765576|emb|Z78451.1|PHZ78451', 'gi|2765610|emb|Z78485.1|PHZ78485', 'gi|2765638|emb|Z78513.1|PBZ78513', 'gi|2765609|emb|Z78484.1|PCZ78484', 'gi|2765594|emb|Z78469.1|PHZ78469', 'gi|2765618|emb|Z78493.1|PGZ78493', 'gi|2765567|emb|Z78442.1|PBZ78442', 'gi|2765622|emb|Z78497.1|PDZ78497', 'gi|2765590|emb|Z78465.1|PRZ78465', 'gi|2765608|emb|Z78483.1|PVZ78483', 'gi|2765605|emb|Z78480.1|PGZ78480', 'gi|2765651|emb|Z78526.1|CGZ78526', 'gi|2765564|emb|Z78439.1|PBZ78439', 'gi|2765593|emb|Z78468.1|PAZ78468', 'gi|2765589|emb|Z78464.1|PGZ78464', 'gi|2765640|emb|Z78515.1|MXZ78515', 'gi|2765604|emb|Z78479.1|PPZ78479', 'gi|2765637|emb|Z78512.1|PWZ78512', 'gi|2765573|emb|Z78448.1|PAZ78448', 'gi|2765636|emb|Z78511.1|PEZ78511', 'gi|2765586|emb|Z78461.1|PWZ78461', 'gi|2765621|emb|Z78496.1|PAZ78496', 'gi|2765623|emb|Z78498.1|PMZ78498', 'gi|2765658|emb|Z78533.1|CIZ78533', 'gi|2765596|emb|Z78471.1|PDZ78471', 'gi|2765627|emb|Z78502.1|PBZ78502', 'gi|2765639|emb|Z78514.1|PSZ78514', 'gi|2765648|emb|Z78523.1|CHZ78523', 'gi|2765625|emb|Z78500.1|PWZ78500', 'gi|2765595|emb|Z78470.1|PPZ78470', 'gi|2765572|emb|Z78447.1|PVZ78447', 'gi|2765614|emb|Z78489.1|PDZ78489', 'gi|2765583|emb|Z78458.1|PHZ78458', 'gi|2765628|emb|Z78503.1|PCZ78503', 'gi|2765652|emb|Z78527.1|CYZ78527', 'gi|2765641|emb|Z78516.1|CPZ78516', 'gi|2765607|emb|Z78482.1|PEZ78482', 'gi|2765617|emb|Z78492.1|PBZ78492', 'gi|2765624|emb|Z78499.1|PMZ78499', 'gi|2765578|emb|Z78453.1|PSZ78453', 'gi|2765592|emb|Z78467.1|PSZ78467', 'gi|2765570|emb|Z78445.1|PUZ78445', 'gi|2765587|emb|Z78462.1|PSZ78462', 'gi|2765629|emb|Z78504.1|PKZ78504', 'gi|2765597|emb|Z78472.1|PLZ78472', 'gi|2765565|emb|Z78440.1|PPZ78440', 'gi|2765602|emb|Z78477.1|PVZ78477', 'gi|2765598|emb|Z78473.1|PSZ78473', 'gi|2765657|emb|Z78532.1|CCZ78532', 'gi|2765635|emb|Z78510.1|PCZ78510', 'gi|2765643|emb|Z78518.1|CRZ78518', 'gi|2765569|emb|Z78444.1|PAZ78444', 'gi|2765574|emb|Z78449.1|PMZ78449', 'gi|2765642|emb|Z78517.1|CFZ78517', 'gi|2765654|emb|Z78529.1|CLZ78529', 'gi|2765603|emb|Z78478.1|PVZ78478', 'gi|2765647|emb|Z78522.1|CMZ78522', 'gi|2765585|emb|Z78460.1|PCZ78460', 'gi|2765591|emb|Z78466.1|PPZ78466', 'gi|2765580|emb|Z78455.1|PJZ78455', 'gi|2765613|emb|Z78488.1|PTZ78488', 'gi|2765646|emb|Z78521.1|CCZ78521', 'gi|2765649|emb|Z78524.1|CFZ78524', 'gi|2765606|emb|Z78481.1|PIZ78481', 'gi|2765619|emb|Z78494.1|PNZ78494', 'gi|2765599|emb|Z78474.1|PKZ78474']

Specifying the dictionary keys

Suppose you want to use the same keys as before? You’ll need to write a tiny function to map from the FASTA identifier (as a string) to the key you want:

In [33]:
def get_acc(identifier):
    """"Given a SeqRecord identifier string, return the accession number as a string.

    e.g. "gi|2765613|emb|Z78488.1|PTZ78488" -> "Z78488.1"
    """
    parts = identifier.split("|")
    assert len(parts) == 5 and parts[0] == "gi" and parts[2] == "emb"
    return parts[3]

Then we can give this function to the Bio.SeqIO.index() function to use in building the dictionary:

In [ ]:
orchid_dict = SeqIO.index("data/ls_orchid.fasta", "fasta", key_function=get_acc)
print(list(orchid_dict.keys()))
['Z78459.1', 'Z78496.1', 'Z78501.1', 'Z78443.1', 'Z78514.1', 'Z78452.1', 'Z78442.1', 'Z78439.1', 'Z78526.1', 'Z78527.1', 'Z78532.1', 'Z78471.1', 'Z78444.1', 'Z78511.1', 'Z78489.1', 'Z78502.1', 'Z78480.1', 'Z78453.1', 'Z78505.1', 'Z78475.1', 'Z78445.1', 'Z78507.1', 'Z78515.1', 'Z78446.1', 'Z78462.1', 'Z78461.1', 'Z78476.1', 'Z78523.1', 'Z78521.1', 'Z78440.1', 'Z78483.1', 'Z78449.1', 'Z78469.1', 'Z78458.1', 'Z78506.1', 'Z78470.1', 'Z78510.1', 'Z78477.1', 'Z78499.1', 'Z78519.1', 'Z78495.1', 'Z78474.1', 'Z78491.1', 'Z78487.1', 'Z78479.1', 'Z78450.1', 'Z78513.1', 'Z78492.1', 'Z78456.1', 'Z78467.1', 'Z78504.1', 'Z78531.1', 'Z78481.1', 'Z78530.1', 'Z78447.1', 'Z78460.1', 'Z78503.1', 'Z78500.1', 'Z78464.1', 'Z78498.1', 'Z78472.1', 'Z78454.1', 'Z78465.1', 'Z78473.1', 'Z78482.1', 'Z78490.1', 'Z78497.1', 'Z78508.1', 'Z78463.1', 'Z78466.1', 'Z78524.1', 'Z78518.1', 'Z78520.1', 'Z78509.1', 'Z78455.1', 'Z78516.1', 'Z78451.1', 'Z78485.1', 'Z78486.1', 'Z78512.1', 'Z78468.1', 'Z78488.1', 'Z78457.1', 'Z78529.1', 'Z78494.1', 'Z78522.1', 'Z78525.1', 'Z78484.1', 'Z78533.1', 'Z78441.1', 'Z78448.1', 'Z78493.1', 'Z78478.1', 'Z78517.1']

Getting the raw data for a record

The dictionary-like object from Bio.SeqIO.index() gives you each entry as a SeqRecord object. However, it is sometimes useful to be able to get the original raw data straight from the file. For this use the get_raw() method which takes a single argument (the record identifier) and returns a string (extracted from the file without modification).

A motivating example is extracting a subset of a records from a large file where either Bio.SeqIO.write() does not (yet) support the output file format (e.g. the plain text SwissProt file format) or where you need to preserve the text exactly (e.g. GenBank or EMBL output from Biopython does not yet preserve every last bit of annotation).

Let’s suppose you have download the whole of UniProt in the plain text SwissPort file format from their FTP site (ftp://ftp.uniprot.org/pub/databases/uniprot/current_release/knowledgebase/complete/uniprot_sprot.dat.gz - Careful big download) and uncompressed it as the file uniprot_sprot.dat, and you want to extract just a few records from it:

[We should find a smaller example here, Tiago]

In [ ]:
#Use this to download the file
!wget -c ftp://ftp.uniprot.org/pub/databases/uniprot/current_release/knowledgebase/complete/uniprot_sprot.dat.gz -O data/uniprot_sprot.dat.gz
!gzip -d data/uniprot_sprot.dat.gz
--2015-12-25 07:24:10--  ftp://ftp.uniprot.org/pub/databases/uniprot/current_release/knowledgebase/complete/uniprot_sprot.dat.gz
           => ‘data/uniprot_sprot.dat.gz’
Resolving ftp.uniprot.org (ftp.uniprot.org)... 141.161.180.197
Connecting to ftp.uniprot.org (ftp.uniprot.org)|141.161.180.197|:21...
In [ ]:
uniprot = SeqIO.index("data/uniprot_sprot.dat", "swiss")
handle = open("data/selected.dat", "w")
for acc in ["P33487", "P19801", "P13689", "Q8JZQ5", "Q9TRC7"]:
    handle.write(uniprot.get_raw(acc).decode("utf-8"))
handle.close()

There is a longer example in sorting section using the SeqIO.index() function to sort a large sequence file (without loading everything into memory at once).

Sequence files as Dictionaries – Database indexed files

Biopython 1.57 introduced an alternative, Bio.SeqIO.index_db(), which can work on even extremely large files since it stores the record information as a file on disk (using an SQLite3 database) rather than in memory. Also, you can index multiple files together (providing all the record identifiers are unique).

The Bio.SeqIO.index() function takes three required arguments:

  • Index filename, we suggest using something ending .idx. This index file is actually an SQLite3 database.
  • List of sequence filenames to index (or a single filename)
  • File format (lower case string as used in the rest of the SeqIO module).

As an example, consider the GenBank flat file releases from the NCBI FTP site, ftp://ftp.ncbi.nih.gov/genbank/, which are gzip compressed GenBank files. As of GenBank release 182, there are 16 files making up the viral sequences, gbvrl1.seq, ..., gbvrl16.seq, containing in total almost one million records. You can index them like this:

In [ ]:
#this will download the files - Currently there are more than 16, but we will do only 4
import os
for i in range(1, 5):
    os.system('wget ftp://ftp.ncbi.nih.gov/genbank/gbvrl%i.seq.gz -O data/gbvrl%i.seq.gz' % (i, i))
    os.system('gzip -d data/gbvrl%i.seq.gz' % i)
In [ ]:
files = ["data/gbvrl%i.seq" % i for i in range(1, 5)]
gb_vrl = SeqIO.index_db("data/gbvrl.idx", files, "genbank")
print("%i sequences indexed" % len(gb_vrl))

That takes about two minutes to run on my machine. If you rerun it then the index file (here gbvrl.idx) is reloaded in under a second. You can use the index as a read only Python dictionary - without having to worry about which file the sequence comes from, e.g.

In [ ]:
print(next(gb_vrl.keys()))
In [ ]:
print(gb_vrl["AB000048.1"].description)

Getting the raw data for a record

Just as with the Bio.SeqIO.index() function, the dictionary like object also lets you get at the raw text of each record:

In [ ]:
print(gb_vrl.get_raw("AB000048.1"))

Indexing compressed files

Very often when you are indexing a sequence file it can be quite large - so you may want to compress it on disk. Unfortunately efficient random access is difficult with the more common file formats like gzip and bzip2. In this setting, BGZF (Blocked GNU Zip Format) can be very helpful. This is a variant of gzip (and can be decompressed using standard gzip tools) popularised by the BAM file format, http://samtools.sourceforge.net/, and http://samtools.sourceforge.net/tabix.shtml.

To create a BGZF compressed file you can use the command line tool bgzip which comes with samtools. In our examples we use a filename extension *.bgz, so they can be distinguished from normal gzipped files (named *.gz). You can also use the Bio.bgzf module to read and write BGZF files from within Python.

The Bio.SeqIO.index() and Bio.SeqIO.index_db() can both be used with BGZF compressed files. For example, if you started with an uncompressed GenBank file:

In [ ]:
orchid_dict = SeqIO.index("data/ls_orchid.gbk", "genbank")
print(len(orchid_dict))

You could compress this (while keeping the original file) at the command line using the following command:

$ bgzip -c ls_orchid.gbk > ls_orchid.gbk.bgz

You can use the compressed file in exactly the same way:

In [ ]:
orchid_dict = SeqIO.index("data/ls_orchid.gbk.bgz", "genbank")
print(len(orchid_dict))

or

In [ ]:
orchid_dict = SeqIO.index_db("data/ls_orchid.gbk.bgz.idx", "data/ls_orchid.gbk.bgz", "genbank")
print(len(orchid_dict))

The SeqIO indexing automatically detects the BGZF compression. Note that you can’t use the same index file for the uncompressed and compressed files.

Discussion

So, which of these methods should you use and why? It depends on what you are trying to do (and how much data you are dealing with). However, in general picking Bio.SeqIO.index() is a good starting point. If you are dealing with millions of records, multiple files, or repeated analyses, then look at Bio.SeqIO.index_db().

Reasons to choose Bio.SeqIO.to_dict() over either Bio.SeqIO.index() or Bio.SeqIO.index_db() boil down to a need for flexibility despite its high memory needs. The advantage of storing the SeqRecord objects in memory is they can be changed, added to, or removed at will. In addition to the downside of high memory consumption, indexing can also take longer because all the records must be fully parsed.

Both Bio.SeqIO.index() and Bio.SeqIO.index_db() only parse records on demand. When indexing, they scan the file once looking for the start of each record and do as little work as possible to extract the identifier.

Reasons to choose Bio.SeqIO.index() over Bio.SeqIO.index_db() include:

  • Faster to build the index (more noticeable in simple file formats)
  • Slightly faster access as SeqRecord objects (but the difference is only really noticeable for simple to parse file formats).
  • Can use any immutable Python object as the dictionary keys (e.g. a tuple of strings, or a frozen set) not just strings.
  • Don’t need to worry about the index database being out of date if the sequence file being indexed has changed.

Reasons to choose Bio.SeqIO.index_db() over Bio.SeqIO.index() include:

  • Not memory limited - this is already important with files from second generation sequencing where 10s of millions of sequences are common, and using Bio.SeqIO.index() can require more than 4GB of RAM and therefore a 64bit version of Python.
  • Because the index is kept on disk, it can be reused. Although building the index database file takes longer, if you have a script which will be rerun on the same datafiles in future, this could save time in the long run.
  • Indexing multiple files together
  • The get_raw() method can be much faster, since for most file formats the length of each record is stored as well as its offset.

Writing sequence files

We’ve talked about using Bio.SeqIO.parse() for sequence input (reading files), and now we’ll look at Bio.SeqIO.write() which is for sequence output (writing files). This is a function taking three arguments: some SeqRecord objects, a handle or filename to write to, and a sequence format.

Here is an example, where we start by creating a few SeqRecord objects the hard way (by hand, rather than by loading them from a file):

In [ ]:
rec1 = SeqRecord(Seq("MMYQQGCFAGGTVLRLAKDLAENNRGARVLVVCSEITAVTFRGPSETHLDSMVGQALFGD" \
                    +"GAGAVIVGSDPDLSVERPLYELVWTGATLLPDSEGAIDGHLREVGLTFHLLKDVPGLISK" \
                    +"NIEKSLKEAFTPLGISDWNSTFWIAHPGGPAILDQVEAKLGLKEEKMRATREVLSEYGNM" \
                    +"SSAC", generic_protein),
                 id="gi|14150838|gb|AAK54648.1|AF376133_1",
                 description="chalcone synthase [Cucumis sativus]")

rec2 = SeqRecord(Seq("YPDYYFRITNREHKAELKEKFQRMCDKSMIKKRYMYLTEEILKENPSMCEYMAPSLDARQ" \
                    +"DMVVVEIPKLGKEAAVKAIKEWGQ", generic_protein),
                 id="gi|13919613|gb|AAK33142.1|",
                 description="chalcone synthase [Fragaria vesca subsp. bracteata]")

rec3 = SeqRecord(Seq("MVTVEEFRRAQCAEGPATVMAIGTATPSNCVDQSTYPDYYFRITNSEHKVELKEKFKRMC" \
                    +"EKSMIKKRYMHLTEEILKENPNICAYMAPSLDARQDIVVVEVPKLGKEAAQKAIKEWGQP" \
                    +"KSKITHLVFCTTSGVDMPGCDYQLTKLLGLRPSVKRFMMYQQGCFAGGTVLRMAKDLAEN" \
                    +"NKGARVLVVCSEITAVTFRGPNDTHLDSLVGQALFGDGAAAVIIGSDPIPEVERPLFELV" \
                    +"SAAQTLLPDSEGAIDGHLREVGLTFHLLKDVPGLISKNIEKSLVEAFQPLGISDWNSLFW" \
                    +"IAHPGGPAILDQVELKLGLKQEKLKATRKVLSNYGNMSSACVLFILDEMRKASAKEGLGT" \
                    +"TGEGLEWGVLFGFGPGLTVETVVLHSVAT", generic_protein),
                 id="gi|13925890|gb|AAK49457.1|",
                 description="chalcone synthase [Nicotiana tabacum]")

my_records = [rec1, rec2, rec3]

Now we have a list of SeqRecord objects, we’ll write them to a FASTA format file:

In [ ]:
SeqIO.write(my_records, "data/my_example.faa", "fasta")

And if you open this file in your favourite text editor it should look like this:

>gi|14150838|gb|AAK54648.1|AF376133_1 chalcone synthase [Cucumis sativus]
MMYQQGCFAGGTVLRLAKDLAENNRGARVLVVCSEITAVTFRGPSETHLDSMVGQALFGD
GAGAVIVGSDPDLSVERPLYELVWTGATLLPDSEGAIDGHLREVGLTFHLLKDVPGLISK
NIEKSLKEAFTPLGISDWNSTFWIAHPGGPAILDQVEAKLGLKEEKMRATREVLSEYGNM
SSAC
>gi|13919613|gb|AAK33142.1| chalcone synthase [Fragaria vesca subsp. bracteata]
YPDYYFRITNREHKAELKEKFQRMCDKSMIKKRYMYLTEEILKENPSMCEYMAPSLDARQ
DMVVVEIPKLGKEAAVKAIKEWGQ
>gi|13925890|gb|AAK49457.1| chalcone synthase [Nicotiana tabacum]
MVTVEEFRRAQCAEGPATVMAIGTATPSNCVDQSTYPDYYFRITNSEHKVELKEKFKRMC
EKSMIKKRYMHLTEEILKENPNICAYMAPSLDARQDIVVVEVPKLGKEAAQKAIKEWGQP
KSKITHLVFCTTSGVDMPGCDYQLTKLLGLRPSVKRFMMYQQGCFAGGTVLRMAKDLAEN
NKGARVLVVCSEITAVTFRGPNDTHLDSLVGQALFGDGAAAVIIGSDPIPEVERPLFELV
SAAQTLLPDSEGAIDGHLREVGLTFHLLKDVPGLISKNIEKSLVEAFQPLGISDWNSLFW
IAHPGGPAILDQVELKLGLKQEKLKATRKVLSNYGNMSSACVLFILDEMRKASAKEGLGT
TGEGLEWGVLFGFGPGLTVETVVLHSVAT

Suppose you wanted to know how many records the Bio.SeqIO.write() function wrote to the handle? If your records were in a list you could just use len(my_records), however you can’t do that when your records come from a generator/iterator. TheBio.SeqIO.write() function returns the number of SeqRecord objects written to the file.

Note - If you tell the Bio.SeqIO.write() function to write to a file that already exists, the old file will be overwritten without any warning.

Round trips

Some people like their parsers to be ‘round-tripable’, meaning if you read in a file and write it back out again it is unchanged. This requires that the parser must extract enough information to reproduce the original file exactly. Bio.SeqIO does not aim to do this.

As a trivial example, any line wrapping of the sequence data in FASTA files is allowed. An identical SeqRecord would be given from parsing the following two examples which differ only in their line breaks:

>YAL068C-7235.2170 Putative promoter sequence
TACGAGAATAATTTCTCATCATCCAGCTTTAACACAAAATTCGCACAGTTTTCGTTAAGA
GAACTTAACATTTTCTTATGACGTAAATGAAGTTTATATATAAATTTCCTTTTTATTGGA

>YAL068C-7235.2170 Putative promoter sequence
TACGAGAATAATTTCTCATCATCCAGCTTTAACACAAAATTCGCA
CAGTTTTCGTTAAGAGAACTTAACATTTTCTTATGACGTAAATGA
AGTTTATATATAAATTTCCTTTTTATTGGA

To make a round-tripable FASTA parser you would need to keep track of where the sequence line breaks occurred, and this extra information is usually pointless. Instead Biopython uses a default line wrapping of \(60\) characters on output. The same problem with white space applies in many other file formats too. Another issue in some cases is that Biopython does not (yet) preserve every last bit of annotation (e.g. GenBank and EMBL).

Occasionally preserving the original layout (with any quirks it may have) is important. See the section about raw access about the get_raw() method of the Bio.SeqIO.index() dictionary-like object for one potential solution.

Converting between sequence file formats

In previous example we used a list of SeqRecord objects as input to the Bio.SeqIO.write() function, but it will also accept a SeqRecord iterator like we get from Bio.SeqIO.parse() - this lets us do file conversion by combining these two functions.

For this example we’ll read in the GenBank format file and write it out in FASTA format:

In [ ]:
records = SeqIO.parse("data/ls_orchid.gbk", "genbank")
count = SeqIO.write(records, "data/my_example.fasta", "fasta")
print("Converted %i records" % count)

Still, that is a little bit complicated. So, because file conversion is such a common task, there is a helper function letting you replace that with just:

In [ ]:
count = SeqIO.convert("data/ls_orchid.gbk", "genbank", "data/my_example.fasta", "fasta")
print("Converted %i records" % count)

The Bio.SeqIO.convert() function will take handles or filenames. Watch out though - if the output file already exists, it will overwrite it! To find out more, see the built in help:

In [ ]:
help(SeqIO.convert)

In principle, just by changing the filenames and the format names, this code could be used to convert between any file formats available in Biopython. However, writing some formats requires information (e.g. quality scores) which other files formats don’t contain. For example, while you can turn a FASTQ file into a FASTA file, you can’t do the reverse.

Finally, as an added incentive for using the Bio.SeqIO.convert() function (on top of the fact your code will be shorter), doing it this way may also be faster! The reason for this is the convert function can take advantage of several file format specific optimisations and tricks.

Converting a file of sequences to their reverse complements

Suppose you had a file of nucleotide sequences, and you wanted to turn it into a file containing their reverse complement sequences. This time a little bit of work is required to transform the SeqRecord objects we get from our input file into something suitable for saving to our output file.

To start with, we’ll use Bio.SeqIO.parse() to load some nucleotide sequences from a file, then print out their reverse complements using the Seq object’s built in .reverse_complement() method:

In [ ]:
for i, record in enumerate(SeqIO.parse("data/ls_orchid.gbk", "genbank")):
    print(record.id)
    print(record.seq.reverse_complement())
    if i == 2:  # 3 is enough
        break

Now, if we want to save these reverse complements to a file, we’ll need to make SeqRecord objects. We can use the SeqRecord object’s built in .reverse_complement() method but we must decide how to name our new records.

This is an excellent place to demonstrate the power of list comprehensions which make a list in memory:

In [ ]:
records = [rec.reverse_complement(id="rc_"+rec.id, description = "reverse complement") \
           for rec in SeqIO.parse("data/ls_orchid.fasta", "fasta")]
len(records)

Now list comprehensions have a nice trick up their sleeves, you can add a conditional statement:

In [ ]:
records = [rec.reverse_complement(id="rc_"+rec.id, description = "reverse complement") \
           for rec in SeqIO.parse("data/ls_orchid.fasta", "fasta") if len(rec)<700]
len(records)

That would create an in memory list of reverse complement records where the sequence length was under 700 base pairs. However, we can do exactly the same with a generator expression - but with the advantage that this does not create a list of all the records in memory at once:

In [ ]:
records = (rec.reverse_complement(id="rc_"+rec.id, description = "reverse complement") \
          for rec in SeqIO.parse("data/ls_orchid.fasta", "fasta") if len(rec)<700)

As a complete example:

In [ ]:
records = (rec.reverse_complement(id="rc_"+rec.id, description = "reverse complement") \
           for rec in SeqIO.parse("data/ls_orchid.fasta", "fasta") if len(rec)<700)
SeqIO.write(records, "data/rev_comp.fasta", "fasta")

Getting your SeqRecord objects as formatted strings

Suppose that you don’t really want to write your records to a file or handle - instead you want a string containing the records in a particular file format. The Bio.SeqIO interface is based on handles, but Python has a useful built in module which provides a string based handle.

For an example of how you might use this, let’s load in a bunch of SeqRecord objects from our orchids GenBank file, and create a string containing the records in FASTA format:

In [ ]:
records = SeqIO.parse("data/ls_orchid.gbk", "genbank")
out_handle = StringIO()
SeqIO.write(records, out_handle, "fasta")
fasta_data = out_handle.getvalue()
print(fasta_data)

This isn’t entirely straightforward the first time you see it! On the bright side, for the special case where you would like a string containing a single record in a particular file format, use the the SeqRecord class’ format() method.

Note that although we don’t encourage it, you can use the format() method to write to a file, for example something like this:

In [ ]:
out_handle = open("data/ls_orchid_long.tab", "w")
for record in SeqIO.parse("data/ls_orchid.gbk", "genbank"):
    if len(record) > 100:
        out_handle.write(record.format("tab"))
out_handle.close()

While this style of code will work for a simple sequential file format like FASTA or the simple tab separated format used here, it will not work for more complex or interlaced file formats. This is why we still recommend using Bio.SeqIO.write(), as in the following example:

In [ ]:
records = (rec for rec in SeqIO.parse("data/ls_orchid.gbk", "genbank") if len(rec) > 100)
SeqIO.write(records, "data/ls_orchid.tab", "tab")

Making a single call to SeqIO.write(...) is also much quicker than multiple calls to the SeqRecord.format(...) method.